hysop.numerics.fft.scipy_fft module¶
FFT iterface for fast Fourier Transforms using scipy fftpack backend (using scipy).
ScipyFFT
ScipyFFTPlan
- class hysop.numerics.fft.scipy_fft.ScipyFFT(backend=None, allocator=None, warn_on_allocation=True, **kwds)[source]¶
Bases:
HostFFTI
Interface to compute local to process FFT-like transforms using the scipy fftpack backend.
- Scipy fftpack backend has some disadvantages:
creation of one or two intermediate temporary buffers at each call
hermitian-complex data layout is different than all other fft implementations
no planning capabilities (scipy.fftpack methods are just wrapped into fake plans)
It has native float and double precision support unlike the numpy fft backend. Planning won’t destroy original inputs.
Initializes the interface and default supported real and complex types.
- dct(a, out=None, type=2, axis=-1, **kwds)[source]¶
Compute the one-dimensional Cosine Transform of specified type.
- Parameters:
a (array_like) – Real input array.
out (array_like) – Real output array of matching input type and shape.
axis (int, optional) – Axis over witch to compute the transform. Defaults to last axis.
- Return type:
(shape, dtype) of the output array determined from the input array.
- dst(a, out=None, type=2, axis=-1, **kwds)[source]¶
Compute the one-dimensional Sine Transform of specified type.
- Parameters:
a (array_like) – Real input array.
out (array_like) – Real output array of matching input type and shape.
axis (int, optional) – Axis over witch to compute the transform. Defaults to last axis.
- Return type:
(shape, dtype) of the output array determined from the input array.
- fft(a, out=None, axis=-1, **kwds)[source]¶
Compute the unscaled one-dimensional complex to complex discrete Fourier Transform.
- Parameters:
a (array_like of np.complex64 or np.complex128) – Complex input array.
out (array_like of np.complex64 or np.complex128) – Complex output array of the same shape and dtype as the input.
axis (int, optional) – Axis over witch to compute the FFT. Defaults to last axis.
- Return type:
(shape, dtype) of the output array determined from the input array.
Notes
N = a.shape[axis] out[0] will contain the sum of the signal (zero-frequency term always real for real inputs).
- If N is even:
out[1:N/2] contains the positive frequency terms. out[N/2] contains the Nyquist frequency (always real for real inputs). out[N/2+1:] contains the negative frequency terms.
- Else if N is odd:
out[1:(N-1)/2] contains the positive frequency terms. out[(N-1)/2:] contains the negative frequency terms.
- idct(a, out=None, type=2, axis=-1, scaling=None, **kwds)[source]¶
Compute the one-dimensional Inverse Cosine Transform of specified type.
- Default scaling is 1/(2*N) for IDCT type (2,3,4) and
1/(2*N-2) for IDCT type 1.
- Parameters:
a (array_like) – Real input array.
out (array_like) – Real output array of matching input type and shape.
axis (int, optional) – Axis over witch to compute the transform. Defaults to last axis.
- Returns:
(shape, dtype, inverse_type, logical_size) of the output array determined from the input
array.
- idst(a, out=None, type=2, axis=-1, scaling=None, **kwds)[source]¶
Compute the one-dimensional Inverse Sine Transform of specified type.
- Default scaling is 1/(2*N) for IDST type (2,3,4) and
1/(2*N+2) for IDST type 1.
- Parameters:
a (array_like) – Real input array.
out (array_like) – Real output array of matching input type and shape.
axis (int, optional) – Axis over witch to compute the transform. Defaults to last axis.
- Returns:
(shape, dtype, inverse_type, logical_size) of the output array determined from the input
array.
- ifft(a, out=None, axis=-1, **kwds)[source]¶
Compute the one-dimensional complex to complex discrete Fourier Transform scaled by 1/N.
- Parameters:
a (array_like of np.complex64 or np.complex128) – Complex input array.
out (array_like of np.complex64 or np.complex128) – Complex output array of the same shape and dtype as the input.
axis (int, optional) – Axis over witch to compute the FFT. Defaults to last axis.
- Return type:
(shape, dtype, logical_size) of the output array determined from the input array.
- irfft(a, out=None, n=None, axis=-1, **kwds)[source]¶
Compute the one-dimensional hermitian complex to real discrete Fourier Transform scaled by 1/N.
- Parameters:
a (array_like of np.complex64 or np.complex128) – Complex input array.
out (array_like of np.float32 or np.float64) –
Real output array of matching real type. out.shape[…] = a.shape[…] Last axis should match forward transform size used:
out.shape[axis] = 2*(a.shape[axis]-1)
out.shape[axis] = 2*(a.shape[axis]-1) + 1
n (int, optional) – Length of the transformed axis of the output. ie: n should be in [2*(a.shape[axis]-1), 2*(a.shape[axis]-1)+1]
axis (int, optional) – Axis over witch to compute the transform. Defaults to last axis.
Notes
To get an odd number of output points, n or out must be specified.
- Returns:
(shape, dtype, logical_size) of the output array determined from the input array,
out and n.
- rfft(a, out=None, axis=-1, **kwds)[source]¶
Compute the unscaled one-dimensional real to hermitian complex discrete Fourier Transform.
- Parameters:
a (array_like of np.float32 or np.float64) – Real input array.
out (array_like of np.complex64 or np.complex128) – Complex output array of matching complex dtype. out.shape[…] = a.shape[…] out.shape[axis] = a.shape[axis]//2 + 1
axis (int, optional) – Axis over witch to compute the transform. Defaults to last axis.
- Return type:
(shape, dtype) of the output array determined from the input array.
Notes
For real inputs there is no information in the negative frequency components that is not already available from the positive frequency component because of the Hermitian symmetry.
N = out.shape[axis] = a.shape[axis]//2 + 1 out[0] will contain the sum of the signal (zero-frequency term, always real). If N is even:
out[1:N/2] contains the positive frequency terms. out[N/2] contains the Nyquist frequency (always real).
- Else if N is odd:
out[1:(N+1)/2] contains the positive frequency terms.
- class hysop.numerics.fft.scipy_fft.ScipyFFTPlan(fn, x, out, axis, scaling=None, **kwds)[source]¶
Bases:
HostFFTPlanI
Wrap a scipy fftpack call (scipy.fftpack does not offer real planning capabilities).
- property input_array¶
Return currently planned input array.
- property output_array¶
Return currently planned output array.